track function
The track function is a core component in Ocelot for simulating particle beam dynamics through a MagneticLattice.
It applies transformations from lattice elements and physics processes (like Space Charge, Wakefields and etc.) to a ParticleArray.
Optionally, it calculates beam envelope (Twiss) parameters at specified intervals using the get_envelope function.
Function Signature
def track(
lattice: MagneticLattice,
p_array: ParticleArray,
navi: Navigator = None,
print_progress: bool = True,
calc_tws: bool = True,
bounds: list = None,
return_df: bool = False,
overwrite_progress: bool = True,
slice: str = None, # Python 'slice' keyword is different
twiss_disp_correction: bool = False
) -> Tuple[Union[List[Twiss], pd.DataFrame], ParticleArray]:
# ... function body ...
# Simplified conceptual flow:
# tws_initial = get_envelope(p_array, ...) if calc_tws else Twiss()
# tws_track = [tws_initial]
# for each step in navigator:
# apply_transfer_maps(p_array, step_maps)
# apply_physics_processes(p_array, step_processes)
# if p_array is empty: break
# if calc_tws:
# current_tws = get_envelope(p_array, ...)
# tws_track.append(current_tws)
# print_progress_if_enabled()
# finalize_physics_processes()
# return tws_track, p_array
pass
Description
The track function simulates the passage of a ParticleArray through a given MagneticLattice. It utilizes a Navigator to determine the sequence of optical element transformations (transfer maps) and physics processes to apply.
- At each step dictated by the navigator, the corresponding transfer maps are applied to the
p_array. - Subsequently, any physics processes scheduled for that step are applied.
- If
calc_twsis true, beam statistics (Twiss parameters, moments, etc.) are calculated from thep_arrayusing theget_envelopefunction after each step. These intermediate Twiss parameters are accumulated and returned. - The tracking continues until the end of the lattice is reached or the particle array becomes empty.
Parameters
-
lattice(MagneticLattice) The magnetic lattice structure through which the particles will be tracked. -
p_array(ParticleArray) The input array of particles to be tracked. This array is modified in place. -
navi(Navigator, optional) The navigator object that orchestrates the tracking steps, including the application of physics processes. IfNone(default), a basicNavigatoris instantiated for the givenlatticewithout any physics processes. Note: If providing a custom navigator, ensure its associated lattice is the same as thelatticeargument. -
print_progress(bool, optional, default:True) IfTrue, prints the tracking progress to the console, indicating the current longitudinal position (z) and the physics processes being applied. -
calc_tws(bool, optional, default:True) IfTrue, Twiss parameters and other beam envelope characteristics are calculated from thep_arrayusingget_envelopeat the beginning and after each significant step in the lattice. The results are collected and returned. IfFalse, this calculation is skipped, and a list of emptyTwissobjects (or an empty DataFrame) is returned for the Twiss part of the output. -
bounds(list, optional, default:None) Passed directly to theget_envelopefunction whencalc_twsisTrue. Defines longitudinal bounds[left_bound, right_bound]in units of (standard deviation ofp_array.tau()) for slicing the particle distribution before calculating Twiss parameters. IfNone, the full bunch is used. -
return_df(bool, optional, default:False) IfTrueandcalc_twsisTrue, the list of calculatedTwissobjects is converted into a pandasDataFramebefore being returned. This can be convenient for analysis and plotting. -
overwrite_progress(bool, optional, default:True) IfTrueandprint_progressisTrue, the progress message will overwrite the previous message on the same line in the console, creating a dynamic progress bar effect. IfFalse, each progress message is printed on a new line. -
slice(str, optional, default:None) Passed directly to theget_envelopefunction whencalc_twsisTrueandboundsare specified. Determines the reference point for longitudinal slicing:None: Uses the mean ofp_array.tau()()."Imax": Uses the longitudinal position of maximum current.
-
twiss_disp_correction(bool, optional, default:False) Passed directly to theauto_dispparameter of theget_envelopefunction whencalc_twsisTrue. IfTrue,get_envelopewill attempt to estimate and subtract the linear dispersion effects from the particle statistics when calculating Twiss parameters. This is useful for obtaining betatron Twiss parameters in dispersive regions.
Returns
A tuple: (tws_results, p_array_final)
-
tws_track(List[Twiss]orpd.DataFrame)- If
calc_twsisTrue: A list ofTwissobjects, each representing the beam envelope parameters at different points along the lattice (including the start). Thescoordinate of eachTwissobject indicates its longitudinal position. Ifreturn_dfisTrue, this will be a pandasDataFramederived from theseTwissobjects. - If
calc_twsisFalse: A list containing initial and final emptyTwissobjects, or an empty DataFrame ifreturn_dfisTrue. (The actual behavior based on your code seems to be: an initial empty Twiss, and subsequent empty Twiss objects for each step). The provided code returns a list containing an initialTwissobject (either calculated or empty) and appends one moreTwissobject pernavi.get_next_step()iteration.
- If
-
p_array_final(ParticleArray) TheParticleArrayafter all tracking steps and physics processes have been applied. This is the samep_arrayobject passed as input, modified in place.
p_arrayThe input p_array object is modified in-place during the execution of the track function. Consequently, the returned p_array_final is the same object as the input p_array, now reflecting its state after tracking.
While returning p_array_final might seem redundant given the in-place modification, this behavior is maintained for consistency with existing scripts and a common functional programming pattern of returning outputs. Users should be aware that any references to the original p_array outside the function will also point to the modified particle data. If the original state of p_array needs to be preserved, create a copy before calling track() (e.g., p_array_copy = p_array.copy()).
Key Operations within track (Conceptual Flow)
-
Initialization:
- A default
Navigatoris created if one isn't provided. - If
calc_twsis enabled, initial Twiss parameters are calculated usingget_envelopewith the specifiedbounds,slice, andtwiss_disp_correctionsettings. This initialTwissobject is stored.
- A default
-
Tracking Loop (Iterating through
navi.get_next_step()): The navigator yields steps, each typically comprising:t_maps: A list of transfer map objects (e.g., from drift, quadrupole).dz: The length of this step.proc_list: A list of physics process objects to be applied during this step.phys_steps: The integration lengths for each process inproc_list.
For each such step from the navigator:
- Apply Transfer Maps: Each transfer map in
t_mapsis applied to thep_array. - Apply Physics Processes: Each physics process in
proc_listis applied to thep_arrayover its correspondingz_step. - Check Particle Count: If
p_array.n(number of particles) becomes zero, tracking stops. - Calculate Twiss (if enabled):
get_envelopeis called again on the (now updated)p_arraywith the samebounds,slice, andtwiss_disp_correctionsettings. Thescoordinate of the resultingTwissobject is updated to reflect the current total path lengthL. This newTwissobject is appended totws_track. - Print Progress (if enabled): Information about the current position
zand applied processes is printed.
-
Finalization:
- After the loop, any
finalize()methods of the physics processes registered with the navigator are called.
- After the loop, any
-
Return Values:
- The accumulated list of
Twissobjects (or DataFrame) and the final state of thep_arrayare returned.
- The accumulated list of
Notes
- The
p_arrayis modified in place throughout the tracking. If you need to preserve the initialp_array, make a copy before callingtrack(e.g.,p_array_copy = p_array.copy()). - The accuracy and types of physical effects included depend on the
Navigatorsetup and thePhysProcobjects added to it. - Parameters
bounds,slice, andtwiss_disp_correctiondirectly influence how the intermediate Twiss parameters are calculated byget_envelope, allowing for flexible analysis (e.g., sliced analysis, dispersion-corrected betatron functions).